The R codes in this script are available for download. This work is licensed under the Creative Commons attribution-noncommercial license. Please share and remix noncommercially, mentioning its origin.
This tutorial aims to help you get started using pomp as a suite of tools for analysis of time series data based on stochastic dynamical systems models. First, we give some conceptual background regarding the class of models—partially observed Markov processes (POMPs)—that pomp handles. We then discuss some preliminaries: installing the package and so on. Next, we show how to simulate a POMP using pomp. We then analyze some data using a few different tools. In so doing, we illustrate some of the package’s capabilities by using its algorithms to fit, compare, and criticize the models using pomp’s algorithms. From time to time, exercises for the reader are given.
As its name implies pomp is focused on partially observed Markov process models. These are also known as state-space models, hidden Markov models, and stochastic dynamical systems models. Such models consist of an unobserved Markov state process, connected to the data via an explicit model of the observation process. We refer to the former as the latent process model (or process model for short) and the latter as the measurement model.
Mathematically, each model is a probability distribution. Let \(Y_n\) denote the measurement process and \(X_n\) the latent state process, then by definition, the process model is determined by the density \(f_{X_n|X_{n-1}}\) and the initial density \(f_{X_0}\). The measurement process is determined by the density \(f_{Y_n|X_n}\). These two submodels determine the full POMP model, i.e., the joint distribution \(f_{X_{0:N},Y_{1:N}}\). If we have a sequence of measurements, \(y^*_{n}\), made at times \(t_n\), \(n=1,\dots,N\), then we think of these data, collectively, as a single realization of the \(Y\) process.
The following schematic shows shows causal relations among the process model, the measurement model, and the data. From the statistical point of view, the key perspective is that the model is, at least hypothetically, the process that generated the data.
Structure of a POMP. Arrows show the direction of causality. The closed loop from the state process to itself indicates the dynamic nature of this Markovian process. Information flow runs in the opposite direction.
Mathematically, a POMP is characterized by two conditions.
These conditions can be represented schematically by the following diagram, which indicates the direct dependencies among model variables.
Conditional independence graph of a POMP. The latent state \(X_n\) at time \(t_n\) is conditionally independent of its history given \(X_{n-1}\). The observation \(Y_n\) is conditionally independent of all other variables given \(X_n\). The underlying time scale can be taken to be either discrete or continuous and the observation times need not be equally spaced.
To implement a POMP model in pomp, we have to specify the measurement and process distributions. Note however, that, for each of the process and the measurement models there are two distinct operations that we might desire to perform:
Following the R convention, we refer to the simulation of \(f_{X_n|X_{n-1}}\) as the rprocess component of the POMP model and the evaluation of \(f_{X_n|X_{n-1}}\) as the dprocess component. Similarly, the simulation of \(f_{Y_n|X_n}\) is the rmeasure component while the evaluation of \(f_{Y_n|X_n}\) is the dmeasure component. Finally, we’ll call a simulator of \(f_{X_0}\) the rinit component. Collectively, we’ll refer to these, and other, similarly basic elements of the model, as the model’s basic components.
A method that makes no use of the dprocess component is said to be “plug-and-play” or to have the “plug-and-play property”. At present, pomp is focused on such methods, so there is no reason to discuss the dprocess component further in this document. In the following, we will illustrate and explain how one specifies the rprocess, rmeasure, and dmeasure components of a model in pomp. We will illustrate this using some simple examples.
To get started, we must install pomp, if it is not already installed. This package cannot yet be downloaded from CRAN (though that will change in the near future). However, the latest version is always available at the package homepage on Github. See the package website for installation instructions.
In this document, we will ultimately learn to implement POMP models using the package’s “C snippet” facility. This allows the user to write model components using snippets of C code, which is then compiled and linked into a running R session. This typically affords a manyfold increase in computation time. It is possible to avoid C snippets entirely by writing model components as R functions, and we will begin by doing so, but the speed-ups afforded by C snippets are typically too good to pass up. To use C snippets, you must be able to compile C codes. Compilers are not by default installed on Windows or Mac systems, so users of such systems must do a bit more work to make use of pomp’s facilities. The installation instructions on the package website give details.
Having dispensed with the preliminaries, we now explore some of the functionality provided by pomp. To assist the reader in following this exploration, the R codes for this document are available.
Let us see how to implement a very simple POMP model. In particular, let’s begin by focusing on the famous Ricker model (Ricker 1954), which posits a nonlinear relationship between the size, \(N(t)\), of a population in year \(t\) and its size, \(N(t+1)\), one year later: \[N(t+1)=r\,N(t)\,\exp\left(1-\frac{N(t)}{K}\right).\tag{1}\] Here, \(r\) and \(K\) are constant parameters, usually termed the intrinsic growth rate and the carrying capacity, respectively. As written, this is a deterministic model: it does not allow for any variability in the population dynamics. Let’s modify the Ricker equation by assuming that \(r\) is not constant, but instead has random variation from year to year. If we assume that the intrinsic growth rate varies from year to year as a lognormal random variable, we obtain \[N(t+1)=r\,N(t)\,\exp\left(1-\frac{N(t)}{K}+\varepsilon(t)\right),\tag{2}\] where \(\varepsilon(t)\sim\mathrm{Normal}(0,\sigma)\). Note that we’ve introduced a new parameter, \(\sigma\), which quantifies the intensity of the noise in the population dynamics. Ecologically speaking, Eq. 2 is a model with environmental stochasticity.
Typically, it is relatively straightforward to simulate a POMP model. To accomplish this in pomp, as we’ve already discussed, we specify the rprocess component of the model. We’ll also need to choose values for the model parameters, \(r\), \(K\), and \(\sigma\). We’ll also need to make a choice regarding the initial condition, \(N(0)\). The simplest choice is to treat \(N(0)=N_0\) as a parameter.
library(pomp)
simulate(t0=0, times=1:20,
params=c(r=1.2,K=200,sigma=0.1,N_0=50),
rinit=function (N_0, ...) {
c(N=N_0)
},
rprocess=discrete_time(
function (N, r, K, sigma, ...) {
eps <- rnorm(n=1,mean=0,sd=sigma)
c(N=r*N*exp(1-N/K+eps))
},
delta.t=1
)
) -> sim1## Warning: 'rmeasure' unspecified: NAs generated.
Notice that we’ve specified the rinit and rprocess components of the model as functions. These functions take as arguments the relevant variables (whether these are state variables or parameters). Importantly, they return named numeric vectors. Names of variables and parameters are very important in pomp. Notice too that we’ve used the discrete_time function, which encodes the fact that our Ricker model is a discrete-time stochastic process (a Markov chain). The first argument of discrete_time is an R function encoding Eq. 2; the second argument specifies the discrete time-step.
Note also that the parameters are furnished in the form of a named vector, and that we’ve specified both t0 and times. The former is the initial time, \(t_0\), i.e., the time at which the initial conditions apply. Since our initial condition is that \(N(0)=N_0\), our initial time is \(t_0=0\). The times argument specifies the observation times \(t_1,\dots,t_N\).
Finally, note that we received a warning about NA values being generated. We will soon see what this is about.
What sort of an object is sim1? If we print it, we see
## <object of class 'pomp'>
sim1 is evidently an object of class ‘pomp’. We refer to these as ‘pomp’ objects.
To get more insight into the structure of sim1, we can use spy:
pomp provides methods for plotting ‘pomp’ objects. For example,
We can also recast a ‘pomp’ object as a data frame:
| time | N |
|---|---|
| 1 | 148.3456 |
| 2 | 274.8027 |
| 3 | 225.5831 |
| 4 | 209.7361 |
| 5 | 248.0677 |
| 6 | 276.0464 |
| 7 | 196.4977 |
| 8 | 248.0970 |
| 9 | 214.0114 |
| 10 | 207.6063 |
| 11 | 258.5272 |
| 12 | 245.7282 |
| 13 | 180.2086 |
| 14 | 234.7697 |
| 15 | 244.2032 |
| 16 | 268.4648 |
| 17 | 199.6789 |
| 18 | 211.7116 |
| 19 | 254.1603 |
| 20 | 222.1691 |
Casting the ‘pomp’ object as a data frame allows us to use ggplot2 graphics:
Tp return to the warning we got when we ran simulate: it was telling us that simulate could not make a random draw from the measurement process because we had not supplied it with any information about this process. In particular, we had not supplied a measurement-model simulator. Let’s now see how to specify the measurement-model simulator, or rmeasure.
Let’s suppose that non-lethal traps are used to sample the population to determine its size. Each year, some number of traps are set out and \(Y_t\) is the number of animals captured. We might posit \[Y_t\;\sim\;\mathrm{Poisson}(b\,N(t)),\] where the parameter \(b\) is proportional to sampling effort, e.g., the number of traps. This is a measurement model, and we can implement a simulator for it by specifying another function:
simulate(t0=0, times=1:20,
params=c(r=1.2,K=200,sigma=0.1,N_0=50,b=0.05),
rinit=function (N_0, ...) {
c(N=N_0)
},
rprocess=discrete_time(
function (N, r, K, sigma, ...) {
eps <- rnorm(n=1,mean=0,sd=sigma)
c(N=r*N*exp(1-N/K+eps))
},
delta.t=1
),
rmeasure=function (N, b, ...) {
c(Y=rpois(n=1,lambda=b*N))
}
) -> sim2Note that, again, the rmeasure function need take only the necessary arguments (and ...) and must return a named numeric vector.
Now, in the preceding chunk of code where we construct sim2, there was some redundancy with our earlier construction of sim1. In particular, we specified the same values of t0, times, rinit, and rprocess as before. Since these specifications were stored in sim1, however, we could have simply added in just the new pieces. For example,
simulate(
sim1,
params=c(r=1.2,K=200,sigma=0.1,N_0=50,b=0.05),
rmeasure=function (N, b, ...) {
c(Y=rpois(n=1,lambda=b*N))
}
) -> sim2As before, we can examine our handiwork:
| time | Y | N |
|---|---|---|
| 1 | 11 | 124.0831 |
| 2 | 13 | 200.2960 |
| 3 | 7 | 267.4025 |
| 4 | 8 | 190.2190 |
| 5 | 12 | 195.8285 |
| 6 | 15 | 293.0790 |
| 7 | 14 | 216.3693 |
| 8 | 10 | 217.2981 |
| 9 | 12 | 269.7404 |
| 10 | 16 | 223.0512 |
| 11 | 11 | 200.4050 |
| 12 | 18 | 255.2297 |
| 13 | 12 | 231.3550 |
| 14 | 12 | 197.5631 |
| 15 | 14 | 263.8953 |
| 16 | 14 | 277.9194 |
| 17 | 15 | 197.6269 |
| 18 | 14 | 292.5811 |
| 19 | 13 | 250.8353 |
| 20 | 12 | 206.4404 |
ggplot(data=gather(
as(sim2,"data.frame"),
variable,value,-time),
aes(x=time,y=value,color=variable))+
geom_line()Notice that now our simulate call produces samples from both the \(N\) and the \(Y\) distributions. simulate will try to sample from the joint distribution of latent states and observables, but will sample from just the latent state process if the rmeasure component is undefined.
We can run multiple simulations using the same parameters:
simulate(sim2,nsim=20) -> sims
ggplot(data=gather(
as.data.frame(sims),
variable,value,Y,N
),
aes(x=time,y=value,color=variable,
group=interaction(.id,variable)))+
geom_line()+
facet_grid(variable~.,scales="free_y")+
labs(y="",color="")We can also run a single simulation from multiple parameters:
p <- parmat(coef(sim2),3)
p["sigma",] <- c(0.05,0.25,1)
colnames(p) <- LETTERS[1:3]
simulate(sim2,params=p,format="data.frame") -> sims
sims <- gather(sims,variable,value,Y,N)
ggplot(data=sims,aes(x=time,y=value,color=variable,
group=interaction(.id,variable)))+
geom_line()+
scale_y_log10()+
expand_limits(y=1)+
facet_grid(variable~.id,scales="free_y")+
labs(y="",color="")In the above, coef extracts the parameter vector stored in sim2. parmat makes a matrix with columns that are copies this vector. The second line modifies this matrix so that the columns are a set of points in parameter space that lie along a line with increasing values of the \(\sigma\) parameter.
We can even run multiple simulations at each one of a set of parameters:
simulate(sim2,params=p,
times=seq(0,3),
nsim=500,format="data.frame") -> sims
ggplot(data=separate(sims,.id,c("parset","rep")),
aes(x=N,fill=parset,group=parset,color=parset))+
geom_density(alpha=0.5)+
# geom_histogram(aes(y=..density..),position="dodge")+
facet_grid(time~.,labeller=label_both,scales="free_y")+
lims(x=c(NA,1000))pomp uses certain syntactic conventions. For example, if x is an object and f is a pomp function, then typically f(x,...) -> y yields an object, y, that contains x plus the results of the f operation, together with additional information pertinent to what f did. Here, ... stands for additional arguments to f. Then, if we call f(y), this operation will perform some kind of repeat of the original operation: though the details will vary with f and y, choices made in the original computation (on x) will typically be respected. To change some of these choices, one can do f(y,...), where the ... stands for options of f that correspond to the new choices.
Morever, since y contains so much information, this information will typically be reused when we apply a different function, g, to y, when it makes sense to do so.
As of version 4.1 R provides a new pipe syntax that is especially convenient for us. It allows us to construct pipelines of R commands, which leads (once one is used to it) to code that is easier to read. In brief, a command in ordinary R style, such as
h(g(f(x,...),...),...)
where f, g, h are R functions, x is some R object, and ... represents additional arguments to each of the functions, must be read from the inside out, which is at odds with the sequential nature of the f → g → h computation. An alternative, of course, is to write something like
y <- f(x,...)
z <- g(y,...)
w <- h(z,...)
which multiplies the entities (y, z, w) one must name and keep track of. In contrast, the new pipeline syntax of R allows us to write
x |> f(...) |> g(...) |> h(...)
The object-oriented structure of pomp makes this kind of programming quite natural.
Some experienced R programmers find the pipeline syntax uncomfortable or unnecessary at first, and debugging pipelined code requires a somewhat different approach. Of course, it is not necessary to adopt this style of programming to use pomp, but it is quite natural and, experience shows, quite addictive!
We will use pipelines freely for the remainder of the document.
Important note: If you opt to load the tidyverse (or individual packages therein), be sure to load pomp after loading these packages. There are some name conflicts between the packages that would otherwise cause pomp functions to be masked.
Modify the sim2 object constructed above to change the measurement model. In particular, assume that \[Y_t\;\sim\;\mathrm{NegBin}(b\,N(t),\theta),\] where \(b\) and \(\theta\) are parameters. (By this notation, we understand that \(Y_t\) is negative binomially distributed with mean \(b\,N(t)\) and variance \(b\,N(t)+(b\,N(t))^2/\theta\).) Make and plot some simulations for \(t=20,\dots,50\) with different values of the \(\theta\) parameter. [Hint: use the rnbinom function with the mu, size parameterization.]
We will illustrate some of the pomp data-analysis algorithms by performing a limited analysis on a set of bird abundance data. The data are from annual censuses of a population of Parus major (Great Tit) in Wytham Wood, Oxfordshire. They were retrieved as dataset #10163 from the Global Population Dynamics Database version 2 (NERC Centre for Population Biology, Imperial College, 2010). The original source is McCleery and Perrins (1991). They are provided as part of the package, in the data frame parus. Here we display the data graphically and in a table:
| year | pop |
|---|---|
| 1960 | 148 |
| 1961 | 258 |
| 1962 | 185 |
| 1963 | 170 |
| 1964 | 267 |
| 1965 | 239 |
| 1966 | 196 |
| 1967 | 132 |
| 1968 | 167 |
| 1969 | 186 |
| 1970 | 128 |
| 1971 | 227 |
| 1972 | 174 |
| 1973 | 177 |
| 1974 | 137 |
| 1975 | 172 |
| 1976 | 119 |
| 1977 | 226 |
| 1978 | 166 |
| 1979 | 161 |
| 1980 | 199 |
| 1981 | 306 |
| 1982 | 206 |
| 1983 | 350 |
| 1984 | 214 |
| 1985 | 175 |
| 1986 | 211 |
The basic data-type provided by pomp—the ‘pomp’ object—is designed as a container for a model and data. Let’s construct such an object with these data, together with the Ricker model we’ve already implemented. We do this with a call to the constructor function pomp:
parus |>
pomp(
times="year", t0=1960,
rinit=function (N_0, ...) {
c(N=N_0)
},
rprocess=discrete_time(
function (N, r, K, sigma, ...) {
eps <- rnorm(n=1,mean=0,sd=sigma)
c(N=r*N*exp(1-N/K+eps))
},
delta.t=1
),
rmeasure=function (N, b, ...) {
c(pop=rpois(n=1,lambda=b*N))
}
) -> rickNotice that the measurement model simulator (rmeasure) is not the same as before, reflecting the fact that the data are named pop, not Y as before. Also notice that we specify the time variable by name, thus allowing the data to dictate the observation times. The t0 argument specifies that the stochastic latent state process is initialized in year 1960.
We’ve already shown how to implement a discrete-time model; let’s now see how to implement a continuous-time model. We’ll implement a very simple model of stable but stochastic population dynamics, the logistic, or Verhulst-Pearl, equation with environmental stochasticity. We’ll write this model as a stochastic differential equation (SDE), specifically an Itô diffusion: \[dN = r\,N\,\left(1-\frac{N}{K}\right)\,dt+\sigma\,N\,dW(t),\tag{3}\] where \(N\) is the population size, \(r\) is a fixed rate, the so-called “Malthusian parameter”, \(K\) is the population’s “carrying capacity”, \(\sigma\) describes the intensity of extrinsic stochastic forcing on the system, and \(dW\) is an increment of a standard Wiener process. [Those unfamiliar with Wiener processes and Itô diffusions will not go far wrong thinking of \(dW(t)\), for each time \(t\), as a normal random variable with mean zero and standard deviation \(\sqrt{dt}\).] To implement this model in pomp, we must tell the package how to simulate this model. The easiest way to simulate such an SDE is via the Euler-Maruyama method. In this approximation, we take a large number of very small steps, each of duration \(\Delta t\). At each step, we hold the right-hand side of the above equation constant, compute \(\Delta N\) using that equation, and increment \(N\) accordingly. pomp gives us the euler function to assist us in implementing the Euler-Maruyama method. To use it, we must encode the computations that take a single step. As before, we can do so by writing a function.
vpstep <- function (N, r, K, sigma, delta.t, ...) {
dW <- rnorm(n=1,mean=0,sd=sqrt(delta.t))
c(N = N + r*N*(1-N/K)*delta.t + sigma*N*dW)
}This function computes the value of \(N(t+\Delta t)\) given the value of \(N(t)\), \(\Delta t\), and the parameters.
We fold this with the data into a ‘pomp’ object via a call to the ‘pomp’-object constructor function, pomp. We’ll also include the same measurement model we used before. Because everything but the ‘rprocess’ component is the same as with the Ricker model, to accomplish this, we can simply do
Notice that we’ve specified an Euler-Maruyama step of about 1 day: it will take 365 of these steps to get us from one observation to the next.
As before, we can plot this ‘pomp’ object, recast it as a data frame, and simulate it for any given choice of the parameters.
The following codes produce several simulations for parameters \(r=0.5\), \(K=2000\), \(\sigma=0.1\) and \(b=0.1\), and plot them on the same axes as the data. Notice that the format="data.frame" and include.data=TRUE options facilitate this. Vary the parameters to try to achieve a better fit to the data, as judged purely “by eye”.
vp |>
simulate(
params=c(r=0.5,K=2000,sigma=0.1,b=0.1,N_0=2000),
format="data.frame", include.data=TRUE, nsim=5) |>
mutate(ds=case_when(.id=="data"~"data",TRUE~"simulation")) |>
ggplot(aes(x=year,y=pop,group=.id,color=ds))+
geom_line()+
labs(color="")To this point, we’ve implemented two models by specifying their basic components as R functions. While this has the virtue of transparency, it puts severe constraints on computational performance due to intrinsic limits in the speed with which R codes can be interpreted. If all we want to do is run a few simulations, this is not a problem. As we’ll see, however, in attempting to perform parameter estimation and other inferences, we will need all the speed we can readily get. We achieve potentially massive speed-ups by implementing our basic model components in a language that can be compiled. pomp makes this easy. The key innovation is the C snippet. Let’s see how to code up the Ricker and Verhulst-Pearl models using C snippets. First, we’ll code up the rinit and rmeasure components, which the two models share.
These are about as simple as one can get. As the name suggests, each is just a snippet of C code: not all variables are declared (in fact, in these examples, none are) and the context of the snippet is not specified. In particular, these snippets are not actually complete C functions, so by themselves, they cannot be compiled. They must not actually violate the rules of C syntax however. Among other things, lines must end with a semicolon (;), variable names must respect C restrictions, etc. Although one can enclose essentially arbitrarily complex C code in a C snippet, one can do quite a lot with very simple snippets. If you are new (even very new) to C, don’t worry: once you master a few simple rules, you’ll be able to code up C snippets just as easily, or even more easily, than you code up R functions and you’ll learn to value the resulting speed-up extremely.
Now we’ll code the Ricker and Verhulst-Pearl simulation steps as C snippets.
Csnippet("
double eps = rnorm(0,sigma);
N = r*N*exp(1-N/K+eps);
") -> rickstepC
Csnippet("
double dW = rnorm(0,sqrt(dt));
N += r*N*(1-N/K)*dt+sigma*N*dW;
") -> vpstepCA few observations are in order. First, note that the local variables eps and dW are declared “double”, the standard C data-type for floating point numbers. Observe that dt is a variable that is defined in the context of the vpstepC snippet; when this snippet is executed, dt will be provided by pomp and will be equal to the size of the Euler step actually being taken. Note also that in each of these snippets, the value of N gets over-written by its new value at the next time-step. This is the goal of the C snippets we supply to specify the ‘rprocess’ component of a model. Finally, notice that neither the state variable N nor any of the parameters are declared. The declarations will be handled in a different way, as we’ll see in a moment.
When furnished with one or more C snippets, pomp will provide the necessary declarations and context, compile the resulting C code, dynamically link to it, and use it whenever the corresponding basic model component is needed. We cause all this to happen when we construct an object of class pomp via a call to the constructor function. Let’s do this for the two models now.
parus |>
pomp(
times="year", t0=1960,
rinit=rinit,
rmeasure=rmeas,
rprocess=discrete_time(rickstepC,delta.t=1),
statenames="N",
paramnames=c("r","K","sigma","b","N_0")
) -> rickCparus |>
pomp(
times="year", t0=1960,
rinit=rinit,
rmeasure=rmeas,
rprocess=euler(vpstepC,delta.t=1/365),
statenames="N",
paramnames=c("r","K","sigma","b","N_0")
) -> vpCIn these calls, we use the statenames and paramnames arguments to indicate which of the undeclared variables in the C snippets rickstepC and vpstepC are state variables and which are fixed parameters. Since dW and eps are declared as local variables within the C snippets themselves, we don’t need to mention them here. The rnorm and rpois functions are part of the R API: see the manual on “Writing R Extensions” for a description of these and the other distribution functions provided as part of the R API. A full set of rules for writing pomp C snippets is given in the package help (?Csnippet).
Using the ‘pomp’ objects just constructed, explore model simulations at a variety of different parameters. As before, plot simulations and data on the same axes for comparison purposes.
Write a C snippet implementing the negative binomial measurement model you explored previously. Fold it into the ‘pomp’ objects just constructed. Remember, there is no need to re-specify components you have already specified: by calling pomp on a ‘pomp’ object you can modify some or all of the basic model components. Plot simulations and data on the same axes, as in the immediately preceding exercise.
As mentioned in the introduction, the dmeasure component is the other side of the rmeasure component. The latter simulates the measurement model whereas the former evaluates the measurement model’s probability density function. The following C snippet encodes the dmeasure component.
Here, dpois again comes from the R API. It takes three arguments, the datum (pop), the Poisson parameter (b*N), and give_log. When give_log=0, dpois returns the Poisson likelihood; when give_log=1, dpois returns the log of this likelihood. When this snippet is executed, pomp will provide the value of give_log according to its needs. It is the user’s responsibility to make sure that the correct value is returned for both possible values of give_log. This is one of the most common places where newbies make mistakes!
Write the dmeasure component for your negative binomial model both as an R function and as a C snippet.
We are now in a position to be able to compute the likelihood of the data given any set of parameters for either of our models. For this purpose, we use the particle filter. This powerful algorithm is at the heart of several of pomp’s inference methods. We won’t describe the theory of the particle filter here. The tutorial by Arulampalam, Maskell, Gordon, and Clapp (2002) explains the theory in an accessible way. The pomp Journal of Statistical Software paper gives pseudocode and some examples. The pomp documentation page lists several other tutorial documents that go into more detail.
In pomp, the simplest version of the particle filter is implemented in the function pfilter. Its only required arguments are the ‘pomp’ object and the number of particles, i.e., the Monte Carlo sample size.
rickC |>
pfilter(Np=1000,
params=c(r=1.2,K=2000,sigma=0.3,N_0=1600,b=0.1),
dmeasure=dmeas,
paramnames="b",statenames="N") -> pfrickNotice that, in this call, we specified the dmeasure component using the C snippet we wrote above. What would have happened had we not specified this?
Notice that, because we here introduced a new C snippet, we again had to indicate which of the undeclared variables in dmeas are parameters and which are latent state variables.
What kind of object is pfrick?
## <object of class 'pfilterd_pomp'>
As a ‘pfilterd_pomp’ object, pfrick contains rickC plus a wealth of information regarding the particle filter operation that created it. For example, we can extract the estimated log likelihood at these (arbitrarily chosen) parameters:
## [1] -148.5748
There is also a plot method for ‘pfilterd_pomp’ objects and one for coercing them to data frames.
| year | pop | ess | cond.logLik |
|---|---|---|---|
| 1960 | 148 | 1000.00000 | -3.879800 |
| 1961 | 258 | 276.91478 | -5.329129 |
| 1962 | 185 | 245.07346 | -5.244873 |
| 1963 | 170 | 183.89624 | -5.551510 |
| 1964 | 267 | 245.35009 | -5.443870 |
| 1965 | 239 | 291.78302 | -5.235775 |
| 1966 | 196 | 274.13378 | -5.152608 |
| 1967 | 132 | 66.23966 | -6.545454 |
| 1968 | 167 | 225.27473 | -5.331978 |
| 1969 | 186 | 249.96622 | -5.263931 |
| 1970 | 128 | 56.39030 | -6.692253 |
| 1971 | 227 | 302.62524 | -5.174603 |
| 1972 | 174 | 202.36717 | -5.441043 |
| 1973 | 177 | 211.26047 | -5.369764 |
| 1974 | 137 | 83.00157 | -6.269213 |
| 1975 | 172 | 238.82351 | -5.239950 |
| 1976 | 119 | 45.15378 | -6.853124 |
| 1977 | 226 | 299.48172 | -5.172145 |
| 1978 | 166 | 184.60553 | -5.521961 |
| 1979 | 161 | 146.71626 | -5.714835 |
| 1980 | 199 | 269.30266 | -5.230779 |
| 1981 | 306 | 193.94797 | -5.754881 |
| 1982 | 206 | 313.84250 | -5.057105 |
| 1983 | 350 | 117.75799 | -6.321934 |
| 1984 | 214 | 296.00247 | -5.170972 |
| 1985 | 175 | 203.08438 | -5.449372 |
| 1986 | 211 | 293.20228 | -5.161943 |
Both of these reveal that pfrick contains information about the effective sample size of the particle filter (ess) and the conditional log likelihood, cond.logLik, which in the notation of the introduction, is \[\log f_{Y_n|Y_{1:n-1}}(y_n^*|y_{1:n-1}^*).\]
The particle filter is a Monte Carlo algorithm. Accordingly, it gives us only a noisy estimate of the likelihood. We can reduce this noise by increasing the number of particles, Np, and we can estimate the magnitude of the Monte Carlo error by running a few independent particle filters. For example:
## [1] -148.4803 -148.5014 -148.8655 -148.6632 -149.2521 -148.1796 -148.0358
## [8] -148.5363 -148.4578 -148.5291
## est se
## -148.5019709 0.1004082
In the first line, notice that we did not need to specify Np, despite the fact that there is no default value of this parameter. Indeed, because pfrick is a ‘pfilterd_pomp’ object, it knows the value of Np that was used in its own creation. By default, this same value is used again when it is passed to pfilter. We could, of course, have used a different value simply by specifying Np in this call to pfilter.
The last line of the preceding code chunk computes the log of the mean of the estimated likelihoods and the standard error on this mean via the delta method. Since the particle filter gives an unbiased estimate of the likelihood (not the log likelihood), this operation is sensible, provided the Monte Carlo error is not too large.
Let’s repeat the operation for the Verhulst-Pearl model, again at arbitrary parameters.
vpC |>
pfilter(Np=1000,dmeasure=dmeas,
params=c(r=0.5,K=2000,sigma=0.1,b=0.1,N_0=2000),
paramnames="b",statenames="N") |>
logLik() |>
replicate(n=10) |>
logmeanexp(se=TRUE) -> ll_vp1
ll_vp1Compute the likelihood for the parameters you found in your attempt to estimate parameters “by eye”.
Trajectory matching is the method of estimating the parameters of a deterministic model by fitting the model to data assuming independent measurement errors. Although pomp’s main focus is on stochastic models, it does provide facilities for trajectory matching. pomp makes a conceptual distinction between the stochastic process and a deterministic skeleton, which we can view as a deterministic model related to the stochastic process’ central tendency. We’ll not go into mathematical details here: instead, we’ll illustrate with two examples.
A deterministic skeleton of the stochastic Ricker model is the Ricker map, Eq. 1. We implement this for pomp, again either as an R function or a C snippet and pass it to pomp functions via the skeleton argument. For example:
rickC |>
pomp(
skeleton=map(
Csnippet("DN = r*N*exp(1-N/K);"),
delta.t=1
),
paramnames=c("r","K"), statenames="N"
) -> rickCHere, the left-hand side of Eq. 1 is indicated by the D prefix: in skeleton snippets, we don’t over-write the state variable N. We indicate that the skeleton is a discrete-time map using the map function.
A deterministic skeleton of the Verhulst-Pearl model is the vectorfield (or ordinary differential equation), \[\frac{dN}{dt} = r\,N\,\left(1-\frac{N}{K}\right).\] We fold this into the vpC ‘pomp’ object so:
vpC |>
pomp(
skeleton=vectorfield(Csnippet("DN = r*N*(1-N/K);")),
paramnames=c("r","K"), statenames="N"
) -> vpCSince the skeleton here is a vectorfield, in this C snippet, DN is filled with the value of the time-derivative of N. See the package help (?Csnippet) for a complete set of rules for writing C snippets.
With the deterministic skeleton in place we can generate trajectories of the skeleton using trajectory. For example:
p <- parmat(c(r=0.5,K=2000,sigma=0.1,b=0.1,N_0=2000),10)
p["N_0",] <- seq(10,3000,length=10)
vpC |>
trajectory(params=p,format="data.frame") |>
ggplot(mapping=aes(x=year,y=N,color=.id,group=.id))+
guides(color="none")+
geom_line()+
theme_bw()As with simulate, one can use trajectory to compute multiple trajectories at once, for varying values of the parameters.
In pomp, the function traj_objfun constructs an objective function quantifying the mismatch between model predictions and data. For this purpose, it uses the dmeasure component of the model. This function can be given to any of the large variety of numerical optimizers available in R and R packages. These optimizers search parameter space to find parameters under which the likelihood of the data, given a trajectory of the deterministic skeleton, is maximized.
We’ll demonstrate using the Verhulst-Pearl model.
vpC |>
traj_objfun(
est=c("K","N_0"),
params=c(r=0.5,K=2000,sigma=0.1,b=0.1,N_0=2000),
dmeasure=dmeas, statenames="N", paramnames="b"
) -> ofunThis invocation of traj_objfun creates an objective function, ofun, that can be used to estimate the three parameters \(K\) and \(N_0\). It will hold the other three parameters, \(r\), \(\sigma\), and \(b\), fixed at the values they are given in params.
Notice that, in this code chunk, we had to specify dmeasure once again. Why? What would have happened had we not done so?
What kind of object is ofun?
## <object of class 'traj_match_objfun'>
‘traj_match_objfun’ objects, like the objective functions created by the pomp functions spect_objfun, probe_objfun, and nlf_objfun, is an R function, but it is also more than an R function. It contains not only the vpC ‘pomp’ object, but it additionally saves information each time it is evaluated. Let’s see how such stateful objective functions make it easy to use a wide range of numerical optimization routines.
We can evaluate ofun at any point in the 2-dimensional \(K\)-\(N_0\) space. For example:
## [1] 1161.543
The value returned by ofun is the negative log likelihood, as returned by the model’s dmeasure component. We can estimate the parameters using, for example, the subplex algorithm implemented in the subplex package:
## $par
## [1] 1968.474 1895.245
##
## $value
## [1] 276.2893
##
## $counts
## [1] 318
##
## $convergence
## [1] 0
##
## $message
## [1] "success! tolerance satisfied"
##
## $hessian
## NULL
Note that fit contains, among other things, the estimated parameters (element par) and the minimized value of the negative log likelihood (value).
To make absolutely certain that ofun remembers these estimates, we evaluate it once at fit$par:
## [1] 276.2893
## r K sigma b N_0
## 0.500 1968.474 0.100 0.100 1895.245
## [1] -276.2893
Then we can, for example, extract the fitted trajectory thus:
ofun |>
trajectory(format="data.frame") |>
ggplot(mapping=aes(x=year,y=N,color=.id,group=.id))+
guides(color="none")+
geom_line()+
theme_bw()We can superimpose the model predictions on the data as follows.
ofun |>
trajectory(format="data.frame") |>
mutate(
pop=coef(ofun,"b")*N,
.id="prediction"
) |>
select(-N) |>
rbind(
parus |>
mutate(
.id="data",
pop=as.double(pop)
)
) |>
spread(.id,pop) |>
ggplot(aes(x=year))+
geom_line(aes(y=prediction))+
geom_point(aes(y=data))+
expand_limits(y=0)+
labs(y="pop")Very commonly, model parameters must obey certain constraints. For example, the parameters in the two models we’ve looked at so far are all constrained to be positive. In estimating parameters, however, one frequently wants to employ a numerical optimization method that does not respect constraints. One way of accomodating such unconstrained optimizers is to transform the parameter space so that the constraints disappear. For example, by log-transforming the \(r\) parameter in the Verhulst-Pearl model (Eq. 3), one obtains the superficially different equation \[dN = e^\rho\,N\,\left(1-\frac{N}{K}\right)\,dt+\sigma\,N\,dW(t),\] where \(\rho=\log{r}\). In this model, \(\rho\) can take any (positive or negative) values while \(r\) remains positive.
We incorporate parameter transformations using the partrans argument to many pomp functions, specifying them using the parameter_trans function. In general, the parameter transformations, like other basic model components, can be supplied using R functions or C snippets. If we are merely log-transforming, logit-transforming, or log-barycentric-transforming parameters, however, it is even easier. The following code chunk implements log transformation of the all the parameters of the Verhulst-Pearl and Ricker models.
We could then provide vpr_partrans as needed to any of the various pomp inference methods, via the partrans argument. For example, to estimate parameters for the Verhulst-Pearl model on the transformed scale, we do
vpC |>
traj_objfun(
est=c("b","K","N_0"),
params=c(r=0.5,K=2000,sigma=0.1,b=0.1,N_0=2000),
partrans=parameter_trans(log=c("K","N_0","b")),
dmeasure=dmeas, statenames="N", paramnames=c("b","K","N_0")
) -> ofun2
subplex(log(c(0.1,2000,1500)),fn=ofun2) -> fit
ofun2(fit$par)## [1] 276.2893
## r K sigma b N_0
## 0.50000000 2506.85237865 0.10000000 0.07852373 2413.59563164
Estimate the parameters for the Verhulst-Pearl model assuming negative binomial errors. Do not attempt, at first, to estimate all parameters simultaneously. Focus on estimating \(K\), \(N_0\) and the parameters of the measurement model. It will probably be helpful to make use of parameter transformations to enforce the model constraints.
Estimate some or all of the parameters of the Ricker model using trajectory matching. It is probably a good idea to use parameter transformations.
Let us now turn to the main focus of the pomp package: parameter estimation for fully stochastic models. Iterated filtering is a method for maximizing the likelihood. The method was introduced in its original form by Ionides, Bretó, and King (2006), and was subsequently much improved by Ionides, Nguyen, Atchadé, Stoev, and King (2015). The latter paper rigorously expounds the theory. An addendum to the pomp *Journal of Statistical Software paper provides pseudocode and a simple example. Finally, a tutorial on the theory and practice is linked from the pomp documentation index. Here, we confine ourselves to demonstrating how the IF2 algorithm (Ionides et al. 2015) is applied to the toy examples we have been discussing.
The search for the maximum likelihood parameter estimates (MLE) is inherently a global one: the likelihood surface can (and often enough does) possess multiple local maxima, any of which might represent a plausible explanation of the data in hand. IF2, however, is a local search algorithm: essentially, it starts at a user-provided “guess” and attempts to move “uphill” toward the nearest local maximum it can find. To solve the global problem, therefore, it is advisable to run IF2 multiple times, starting from different “guesses”. By studying the results of such a calculation, one hopes to obtain not only the global MLE but, more importantly, an understanding of the shape of the likelihood surface itself.
To illustrate this, let us attempt to determine at which values of the four parameters \(r\), \(K\), \(\sigma\), and \(N_0\) our Verhulst-Pearl model best explains the Parus major data. Because this is only a toy problem, we will restrict ourselves to the assumption that \(b=1\), despite the fact that it might plausibly be otherwise. In a real scientific study, of course, we would be much less cavalier about fixing parameters arbitrarily in this way.
We begin by constructing a large number of “guesses”. Our design is to distribute these around some multidimensional “box” that, we hope, contains the MLE. Importantly, although our hope is that the box we construct here embraces the MLE, our eventual success does not necessarily depend on it doing so: the IF2 algorithm will be free to search parameter space outside the box. We use the sobol_design function to construct an overdispersed sequence of points within the selected box.
sobol_design(
lower=c(r=0,K=100,sigma=0,N_0=150,b=1),
upper=c(r=5,K=600,sigma=2,N_0=150,b=1),
nseq=100
) -> guesses
plot(guesses,pch=16)In the above call, lower and upper give the bounds of the box and nseq specifies the number of guesses we desire. Now, from each of these starting guesses, we will run the IF2 algorithm (implemented as mif2 in pomp). As its name suggests, this algorithm works by iteratively applying the particle filter: Crucially, IF2 adds an extra ingredient to the particle filter. As the filter is applied, IF2 adds random perturbations to the parameters. This serves several purposes, including smoothing the likelihood surface, combating particle depletion, and allowing IF2 to obtain information about the local structure of the surface. However, this augmented stochastic variability represents an alteration of the model. IF2 therefore gradually reduces the intensity of the random perturbations.
Let us examine a single call to mif2.
vpC |>
mif2(
params=guesses[1,],
Np=1000,
Nmif=20,
dmeasure=dmeas,
partrans=parameter_trans(log=c("r","K","sigma","N_0")),
rw.sd=rw_sd(r=0.02,K=0.02,sigma=0.02,N_0=ivp(0.02)),
cooling.fraction.50=0.5,
paramnames=c("r","K","sigma","N_0","b"),
statenames=c("N")
) -> mf1In this call, we specify the starting point with params: here, we’ve just taken this to be the first of the guesses. As in pfilter, Np sets the number of particles to be used. Nmif fixes the number of IF2 iterations that will be performed. We furnish the necessary dmeasure component via the dmeasure argument, and provide parameter transformations via partrans to enforce positivity of the model parameters. The rw.sd argument specifies the intensity of the perturbations that will be applied at each observation time. The perturbations are normally distributed on the transformed (here, the log) scale, so these represent roughly 2% perturbations per observation. The ivp syntax indicates that N_0 is an initial value parameter and accordingly should be treated slightly differently than the regular parameters. The cooling.fraction.50 argument controls the speed at which the magnitude of the perturbations is drawn down. Finally, the paramnames and statenames arguments are needed since dmeas and the parameter transformations are built using C snippets.
What kind of object is mf1?
## <object of class 'mif2d_pomp'>
pomp provides various methods for ‘mif2d_pomp’ objects. The plot method, for example, produces a diagnostic plot.
The first plot (“filter diagnostics”) shows the results of the last iteration of the particle filter, comparable to the plot produced by calling plot on a ‘pfilterd_pomp’ object. The second plot (“convergence diagnostics”) summarizes the results of the sequence of IF2 iterations. We observe that the log likelihood has increased, concomitant with movement in the \(r\), \(K\), and \(\sigma\) paramters and, to a lesser extent, the \(N_0\) parameter. We observe that, as intended, the \(b\) parameter has remained fixed.
Now, although we have observed the intended improvement in the log likelihood, we should be careful to note that the log likelihood displayed in this plot is the log likelihood of the perturbed model. This model differs from the one we are interested in. To compute the likelihood of our focal model at the parameter returned by mif2, we need to perform a few particle filter operations:
## est se
## -143.9880041 0.1035849
Now we repeat this calculation from each of our 100 guesses. Because these computations are each independent of the others, it is both easy and worthwhile to parallelize them. The foreach package provides some useful tools in this regard.
library(foreach)
foreach (guess=iter(guesses,"row"),
.combine=c, .packages=c("pomp"),
.errorhandling="remove", .inorder=FALSE) %dopar% {
mf1 |> mif2(params=guess)
} -> mifsLet us now walk through the above code. At each of the guesses (foreach (guess=), which are stored in the rows of the guesses data frame (iter(guesses,"row")), we run 20 iterations of the IF2 algorithm (mif2). Notice that the settings of mif2 are not explicitly given here. They are inherited from mf1, which we computed previously. If we had wanted to modify any of them (e.g., Nmif, Np, rw.sd, etc.), we would have simply done so in the call to mif2. Next, we run a few independent particle filters to estimate the log likelihood. Finally, we return the result as a one-row data frame containing estimates of the parameters (coef(mf)), the log likelihood, and its standard error.
We can examine the results of this computation. The following code chunk extracts the traces of the IF2 calculations. These are the movements of the parameters and the likelihood plotted against filter iteration.
mifs |>
traces() |>
melt() |>
filter(variable!="b") |>
ggplot(aes(x=iteration,y=value,group=L1,color=L1))+
geom_line()+
facet_wrap(~variable,scales="free_y")+
guides(color="none")Next, we examine the diagnostics produced in the last iteration of the filter.
mifs |>
as("data.frame") |>
gather(variable,value,-year,-.id) |>
ggplot(aes(x=year,y=value,group=.id,color=.id))+
geom_line()+
facet_wrap(~variable,scales="free_y",ncol=1)+
guides(color="none")As we’ve discussed, to properly evaluate the likelihood at these estimates, we must perform several independent particle filter operations at each one:
foreach (mf=mifs,
.combine=rbind, .packages=c("pomp"),
.errorhandling="remove", .inorder=FALSE) %dopar% {
mf |> pfilter() |> logLik() |> replicate(n=5) |>
logmeanexp(se=TRUE) -> ll
data.frame(as.list(coef(mf)),loglik=ll[1],loglik.se=ll[2])
} -> estimatesA useful way of displaying the results of such multi-start search is the pairwise scatterplot matrix:
estimates |>
full_join(guesses) |>
filter(is.na(loglik) | loglik>max(loglik,na.rm=TRUE)-30) |>
{\(dat)
pairs(~loglik+r+sigma+K+N_0,data=dat,pch=16,
col=if_else(is.na(dat$loglik),"#99999955","#ff0000ff"))
}()The above foreach call, as written, will result in serial, not parallel, execution of the 100 computations. To actually parallelize, one must register a parallel backend. A good first choice is doParallel, which is a separate package that you must explicitly load. The following chunk of code shows how to do this, and prints the number of “workers” that result.
doParallel allows one to exploit multiple cores on a single machine, or to expand a parallel computation across a cluster.
Use doParallel to parallelize the estimation of likelihood for your negative-binomial model.
At this point, there’s no particular reason to suspect that our IF2 searches have arrived at their destination. In general, it is hard to know a priori how much effort will be required to find the MLE. Let’s continue the search, starting with the best points we’ve uncovered so far.
estimates |>
filter(!is.na(loglik)) |>
filter(loglik > max(loglik)-30) |>
select(-loglik,-loglik.se) -> starts
foreach (start=iter(starts,"row"),
.combine=rbind, .packages=c("pomp"),
.errorhandling="remove", .inorder=FALSE) %dopar% {
mf1 |>
mif2(params=start) |>
mif2() -> mf
mf |> pfilter() |> logLik() |> replicate(n=5) |>
logmeanexp(se=TRUE) -> ll
data.frame(as.list(coef(mf)),loglik=ll[1],loglik.se=ll[2])
} -> ests1Note that, in the above, we’ve performed a total of 100 mif2 iterations per starting point. The code above also does the post-mif2 likelihood estimation. It returns just the parameter and likelihood estimates.
We can combine the new estimates with the old ones into a general database:
| r | K | sigma | N_0 | b | loglik | loglik.se | |
|---|---|---|---|---|---|---|---|
| se210 | 3.073142 | 215.2133 | 0.5704475 | 149.9442 | 1 | -141.5969 | 0.1797394 |
| se53 | 3.275665 | 210.0892 | 0.5893403 | 148.5220 | 1 | -141.6427 | 0.1724702 |
| se441 | 3.440774 | 212.4449 | 0.6296984 | 144.9650 | 1 | -141.6554 | 0.2350931 |
| se8 | 1.769714 | 210.7772 | 0.4440766 | 150.4294 | 1 | -141.6664 | 0.0779041 |
| se67 | 1.943607 | 216.3494 | 0.4931726 | 148.8008 | 1 | -141.7087 | 0.1475773 |
| se90 | 5.564673 | 209.7477 | 0.8531969 | 151.0701 | 1 | -141.7094 | 0.2048600 |
estimates |>
filter(loglik>max(loglik,na.rm=TRUE)-4) |>
{\(dat)
pairs(~loglik+r+sigma+K+N_0,data=dat,pch=16)
}()In this plot, we begin to see the emergence of structure in the likelihood surface. In particular, what looks like a ridge of high likelihood is visible in the \(r\)-\(\sigma\) projection. Such structures are very interesting in that they contain clues as to the manner in which the model is fitting the data. They can also pose challenges to efficient estimation, since climbing up to a ridge is harder than traversing it.
We can improve the quality of our estimates and obtain likelihood-ratio-test-based confidence intervals by constructing profile likelihoods. In a likelihood profile, one varies the focal parameter (or parameters) across some range, maximizing the likelihood over the remaining parameters at each value of the focal parameter. The following codes construct a likelihood profile over \(r\).
estimates |>
filter(loglik>max(loglik)-10) |>
select(r,K,sigma,N_0,b) |>
apply(2,range) -> ranges
ranges## r K sigma N_0 b
## [1,] 0.02862099 201.0567 0.3136465 137.9456 1
## [2,] 23.51068046 499.8607 1.4658081 161.2613 1
profile_design(
r=10^seq(
from=log10(ranges[1,1]),
to=log10(ranges[2,1]),
length=20
),
lower=ranges[1,-1],
upper=ranges[2,-1],
nprof=50
) -> starts
dim(starts)## [1] 1000 5
foreach (start=iter(starts,"row"),
.combine=rbind, .packages=c("pomp"),
.errorhandling="remove", .inorder=FALSE) %dopar% {
mf1 |>
mif2(
params=start,
partrans=parameter_trans(log=c("K","sigma","N_0")),
rw.sd=rw_sd(K=0.02,sigma=0.02,N_0=ivp(0.02)),
paramnames=c("K","sigma","N_0","b")
) |>
mif2() -> mf
mf |> pfilter() |> logLik() |> replicate(n=5) |>
logmeanexp(se=TRUE) -> ll
data.frame(as.list(coef(mf)),loglik=ll[1],loglik.se=ll[2])
} -> r_profNotice that we’ve changed the mif2 perturbations (rw.sd): we’ve removed the perturbation on the \(r\) parameter, since we want to hold this parameter fixed.
We add these points to our database:
Next, we plot the likelihood profile. The following plot shows the top two estimates for each value of \(r\) with error bars showing ±2 s.e. and a loess smooth.
r_prof |>
group_by(r) |>
filter(rank(-loglik)<=2) |>
ungroup() |>
ggplot(aes(x=r,y=loglik,
ymin=loglik-2*loglik.se,ymax=loglik+2*loglik.se))+
geom_point()+
geom_errorbar()+
geom_smooth(method="loess",span=0.2)+
scale_x_log10()We see that the 95% likelihood ratio test confidence interval (CI) appears to be one-sided: the CI is roughly \(r>0.75\).
If we plot the profile trace of \(\sigma\), we see that, as we increase \(r\), we have to increase the intensity of the environmental stochasticity to maintain a good fit. Why is this?
r_prof |>
group_by(r) |>
filter(rank(-loglik)<=2) |>
ungroup() |>
ggplot(aes(x=r,y=sigma))+
geom_point()+
geom_smooth(method="loess",span=0.2)+
scale_x_log10()+
labs(y=expression(sigma))Construct a likelihood profile for the \(K\) parameter. Plot the profile and profile traces. Comment on your findings.
r_prof |>
group_by(r) |>
filter(rank(-loglik)<=2) |>
ungroup() |>
select(-loglik,-loglik.se) |>
pivot_longer(-r) |>
group_by(name) |>
summarize(value=predict(loess(value~r),newdata=data.frame(r=2))) |>
ungroup() |>
pivot_wider() |>
bind_cols(r=2) -> theta
bake(file="hindcast1.rds",{
registerDoRNG(174423157)
foreach (i=1:200) %dopar% {
mf1 |> pfilter(params=theta,Np=500,filter.traj=TRUE) -> pf
list(loglik=logLik(pf),traj=filter_traj(pf))
} -> fts
}) -> fts
fts |>
sapply(getElement,"loglik") -> ll
logmeanexp(ll,ess=TRUE)## est ess
## -141.9050 151.6353
fts |>
lapply(getElement,"traj") |>
melt() |>
rename(.id=L1) |>
select(-rep) |>
left_join(
tibble(
.id=as.character(seq_along(ll)),
loglik=ll
),
by=".id"
) |>
mutate(
year=time(mf1)[time]
) |>
group_by(variable,year) |>
summarize(
label=c("lo","med","hi"),
p=c(0.05,0.5,0.95),
q=wquant(value,weights=exp(loglik-max(loglik)),probs=p)
) |>
ungroup() |>
select(-p) |>
pivot_wider(names_from=label,values_from=q) -> quants1
quants1 |>
ggplot()+
geom_line(aes(x=year,y=med),color="darkblue")+
geom_ribbon(aes(x=year,ymin=lo,ymax=hi),color=NA,fill="lightblue",alpha=0.5)+
geom_point(data=parus,aes(x=year,y=pop))+
labs(y="N")pmcmc chainsIf we seek to do full-information Bayesian inference, we can use particle Markov chain Monte Carlo, implemented in pomp in the pmcmc function. The following codes cause parallel pMCMC chains to be run, each beginning at one of the estimates obtained using mif2. Note that we add a prior by furnishing a prior probability density evaluation function (dprior=...). In this case, we assume a product prior: marginally uniform in each of the parameters \(r\), \(\sigma\), \(K\), and \(N_0\).
r_prof |>
group_by(r) |>
filter(loglik==max(loglik)) |>
ungroup() |>
filter(
r > 0.75,
r < 6,
sigma < 2,
K > 100,
K < 600,
N_0 > 100,
N_0 < 600
) |>
select(-loglik,-loglik.se) -> starts
starts| r | K | sigma | N_0 | b |
|---|---|---|---|---|
| 0.978756 | 215.4200 | 0.3669298 | 151.7288 | 1 |
| 1.393392 | 209.7336 | 0.3978699 | 152.4494 | 1 |
| 1.983683 | 204.1178 | 0.4818092 | 150.3479 | 1 |
| 2.824043 | 208.0662 | 0.5682956 | 147.6887 | 1 |
| 4.020410 | 212.5271 | 0.6613447 | 147.5343 | 1 |
| 5.723600 | 212.1773 | 0.7927863 | 153.0553 | 1 |
foreach (start=iter(starts,"row"),.combine=c,
.errorhandling="remove",.inorder=FALSE) %dopar%
{
library(pomp)
mf1 |>
pmcmc(
Nmcmc=2000,Np=200,params=start,
dprior=Csnippet("
lik = dunif(r,0,10,1)+dunif(sigma,0,2,1)+
dunif(K,0,600,1)+dunif(N_0,0,600,1);
lik = (give_log) ? lik : exp(lik);"
),
paramnames=c("K","N_0","r","sigma"),
proposal=mvn_rw_adaptive(
rw.sd=c(r=0.02,sigma=0.02,K=50,N_0=50),
scale.start=100,shape.start=100
)
) -> chain
chain |>
pmcmc(
Nmcmc=20000,
proposal=mvn_rw(covmat(chain))
)
} -> chainsIn the above, we first subset out those mif2 estimates that are consistent with our prior. At each of these, we then perform a number of MCMC moves, using an adaptive multivariate-normal random-walk proposal distribution mvn.rw.adaptive. To complete the inference, we do more MCMC iterations using a multivariate random-walk proposal (mvn.rw) with covariance matrix derived from the preceding computation (covmat).
Now we investigate the chains, to try and determine whether they have converged.
## r sigma K N_0
## 0.8072667 0.8072667 0.8072667 0.8072667
We see that the rejection rate is very good. Let us examine the autocorrelation in the chains.
## loglik log.prior r K sigma N_0 b
## Lag 1 0.92190985 NaN 0.95461111 0.89013855 0.93779316 0.915191693 NaN
## Lag 5 0.68450312 NaN 0.80773770 0.58415560 0.74022358 0.650342606 NaN
## Lag 10 0.49070945 NaN 0.67556882 0.36633158 0.57138932 0.437747486 NaN
## Lag 50 0.09675077 NaN 0.21996975 0.02305981 0.12242624 0.050617555 NaN
## Lag 100 0.01789541 NaN 0.07103979 0.01015425 0.02113539 0.009386808 NaN
## loglik log.prior r K sigma N_0 b
## 4241.592 0.000 2337.883 5962.378 3493.462 5049.668 0.000
## loglik log.prior r K sigma N_0 b
## 1032.086 0.000 1059.283 1060.125 1339.431 1099.495 0.000
The autocorrelation is strong, but drops to small values by lag 100. Accordingly, we thin the chains by a factor of 100. We discard 2000 iterations as a burn-in.
Now let us examine the traces.
traces |>
lapply(as.data.frame) |>
lapply(rownames_to_column,"iter") |>
bind_rows(.id="chain") |>
mutate(iter=as.numeric(iter)) |>
select(chain,iter,loglik,r,sigma,K,N_0) |>
pivot_longer(c(-chain,-iter)) |>
ggplot(aes(x=iter,group=chain,color=chain,y=value))+
guides(color="none")+
labs(x="iteration",y="")+
geom_line(alpha=0.3)+geom_smooth(method="loess",se=FALSE)+
facet_wrap(name~.,scales="free_y",strip.position="left",ncol=2)+
theme(
strip.placement="outside",
strip.background=element_rect(fill=NA,color=NA)
)## Potential scale reduction factors:
##
## Point est. Upper C.I.
## r 1.004 1.018
## sigma 1.002 1.010
## K 0.998 0.999
## N_0 1.000 1.006
##
## Multivariate psrf
##
## 1.01
The trace-plots show good mixing and the Gelman-Rubin statistic confirms that the chains are mixing among themselves: a good indicator of convergence. Thus, insofar as we can tell by these diagnostics, the MCMC iterations appear to have converged to their stationary distribution and we are sampling it well after thinning.
The plots of the marginal posterior densities are shown below.
traces |>
lapply(as.data.frame) |>
lapply(rownames_to_column,"iter") |>
bind_rows(.id="chain") |>
select(chain,iter,loglik,r,sigma,K,N_0) |>
pivot_longer(c(-chain,-iter)) |>
ggplot(aes(x=value))+
geom_density()+
geom_rug()+
labs(x="")+
facet_wrap(~name,scales="free",strip.position="bottom")+
theme(
strip.placement="outside",
strip.background=element_rect(fill=NA,color=NA)
)##
## Iterations = 2000:20000
## Thinning interval = 100
## Number of chains = 6
## Sample size per chain = 181
##
## 1. Empirical mean and standard deviation for each variable,
## plus standard error of the mean:
##
## Mean SD Naive SE Time-series SE
## loglik -142.9970 1.5116 0.045869 0.047658
## log.prior -15.7896 0.0000 0.000000 0.000000
## r 6.0714 2.5371 0.076987 0.087108
## K 213.7084 12.8094 0.388699 0.396234
## sigma 0.8576 0.2353 0.007142 0.007116
## N_0 149.2178 11.9064 0.361298 0.363893
## b 1.0000 0.0000 0.000000 0.000000
##
## 2. Quantiles for each variable:
##
## 2.5% 25% 50% 75% 97.5%
## loglik -146.5687 -143.8087 -142.7488 -141.982 -140.710
## log.prior -15.7896 -15.7896 -15.7896 -15.790 -15.790
## r 1.4522 3.9910 6.3179 8.255 9.835
## K 191.0026 205.3645 212.5041 220.964 242.402
## sigma 0.4231 0.6859 0.8592 1.024 1.318
## N_0 126.7383 141.0953 149.4083 157.208 172.096
## b 1.0000 1.0000 1.0000 1.000 1.000
traces |>
lapply(as.data.frame) |>
lapply(rownames_to_column,"iter") |>
bind_rows(.id="chain") |>
{\(dat)
pairs(~r+sigma+K+N_0,data=dat,pch=16)
}()chains |>
filter_traj() |>
melt() |>
filter(rep > 1000, rep %% 100 == 0) |>
mutate(year=time(mf1)[time]) |>
pivot_wider(names_from=variable) |>
group_by(year) |>
summarize(
label=c("lo","med","hi"),
p=c(0.025,0.5,0.975),
q=wquant(N,probs=p)
) |>
ungroup() |>
select(-p) |>
pivot_wider(names_from=label,values_from=q) -> quants2
quants2 |>
ggplot()+
geom_line(aes(x=year,y=med),color="darkblue")+
geom_ribbon(aes(x=year,ymin=lo,ymax=hi),color=NA,fill="lightblue",alpha=0.5)+
geom_point(data=parus,aes(x=year,y=pop))+
labs(y="N")bind_rows(
with=quants2,
without=quants1,
.id="uncert"
) |>
ggplot()+
geom_line(aes(x=year,y=med,color=uncert))+
geom_ribbon(aes(x=year,ymin=lo,ymax=hi,fill=uncert),color=NA,alpha=0.4)+
geom_point(data=parus,aes(x=year,y=pop))+
labs(y="N")Estimating model parameters by fitting the model to data is typically only a step in the process of trying to understand the processes that generated the data. The next step involves trying to understand how and why the model fits the data the way it does, whether it fits it well, and what scope for improvement there might be. pomp provides a number of tools to facilitate answering these questions through interaction with a fitted model.
Ultimately, since the model is viewed, at least hypothetically, as the process that generated the data, simulation of the fitted model is a central tool we have for model criticism. Let’s plot the data and several simulated realizations of the model process on the same axes.
r_prof |>
filter(loglik==max(loglik)) -> mle
mlepomp <- as(mifs[[1]],"pomp")
coef(mlepomp) <- mle
mlepomp |>
simulate(nsim=8,format="data.frame",include.data=TRUE) |>
ggplot(mapping=aes(x=year,y=pop,group=.id,alpha=(.id=="data")))+
scale_alpha_manual(values=c(`TRUE`=1,`FALSE`=0.2),
labels=c(`FALSE`="simulation",`TRUE`="data"))+
labs(alpha="")+
geom_line()+
theme_bw()The first lines above simply extract the maximum likelihood estimates (mle) from our profle computation. The next pair of lines plug these MLE parameters into a ‘pomp’ object (mlepomp) containing the model and the data. The last set of lines do the simulation and the plotting.
Although it is clear from these plots that the estimated model has more variability and is thus able to explain the data better, it can be hard to read much from spaghetti plots such as this. It’s almost always a good idea to plot the data together with several simulated realizations in order to help assess how similar the two are.
mlepomp |>
simulate(nsim=11,format="data.frame",include.data=TRUE) |>
ggplot(mapping=aes(x=year,y=pop,group=.id,color=(.id=="data")))+
scale_color_manual(values=c(`TRUE`="black",`FALSE`="grey50"),
labels=c(`FALSE`="simulation",`TRUE`="data"))+
labs(color="")+
geom_line()+
facet_wrap(~.id)+
theme_bw()Visual comparison of simulations and data is always a good idea. An indication that the data are not a plausible realization of the model is evidence for lack of fit. In particular, if we have any set of summary statistics, or probes, we can apply them to both simulated and actual data. The probe function facilitates this comparison. Let’s perform this operation using several of the summary statistics provided with pomp: we’ll use the mean, several quantiles, and the autocorrelation at lags 1 and 3.
mlepomp |>
probe(nsim=200,probes=list(
mean=probe_mean("pop"),
q=probe_quantile("pop",probs=c(0.05,0.25,0.5,0.75,0.95)),
probe_acf("pop",lags=c(1,3),type="corr",transform=log)
)) -> vp_probe
vp_probe## <object of class 'probed_pomp'>
For ‘probed_pomp’ object, there are summary and plot methods. There is also an as.data.frame method.
## $coef
## r K sigma N_0 b loglik
## 1.9836835 204.1177581 0.4818092 150.3478525 1.0000000 -141.4331763
## loglik.se
## 0.1645381
##
## $nsim
## [1] 200
##
## $quantiles
## mean q.5% q.25% q.50% q.75% q.95% acf[1] acf[3]
## 0.650 0.595 0.780 0.415 0.500 0.815 0.380 0.840
##
## $pvals
## mean q.5% q.25% q.50% q.75% q.95% acf[1] acf[3]
## 0.6965174 0.8159204 0.4278607 0.8358209 0.9850746 0.3781095 0.7661692 0.3283582
##
## $synth.loglik
## [1] -19.00884
Evidently, summary returns a list with several elements. The quantiles element contains, for each probe, what fraction of the nsim simulations had probe values below the value of the probe applied to the data. The pvals element contains \(P\)-values associated with the two-sided test of the hypothesis that the data were generated by the model.
The plot depicts the multivariate distribution of the probes under the model, with the data-values superimposed. On the diagonal, we see the marginal distributions of the individual probes, represented as histograms, with the vertical line signifying the value of the corresponding probe on the data. Above the diagonal, the scatterplots show the pairwise distributions of probes and the crosshairs, the corresponding data-values. Below the diagonal, the panels contains the pairwise correlations among the simulated probes.
To this point, we’ve seen how to implement POMP models, simulate them, to compute and maximize likelihoods, and perform certain kinds of diagnostic checks. The pomp website contains more documentation, including the full package manual, and a variety of tutorials and short courses. The package itself contains a number of built-in examples and datasets that can be explored.
pomp provides a large toolbox of different inference methods, only a few of which have been explored here. In particular, the package provides other methods for parameter estimation, both in the frequentist and Bayesian modes. See for example (abc, pmcmc, probe.match, spect.match, nlf, enkf, bsmc2). It also provides a variety of tools for model checking (spect, nlf). It is frequently the case that an approach that makes use of more than one approach has advantages over more “purist” approaches: the main goal of the package is to facilitate effective inference by bringing a variety of tools, with complementary strengths and weaknesses, to the user in a common format.
Although the goal of this document has been to introduce the beginner to the package through a display of the pomp toolbox in the context of a rudimentary and incomplete data analysis of a short time series with toy models, it is important to realize that these tools have proven their utility on some extremely challenging problems, including some for which other existing methods are either less efficient or entirely infeasible. The bibliography has links to peer-reviewed publications that have used these methods.
This document was produced with the following software versions:
| R | 4.2.2 |
| pomp | 4.5.4.0 |
| coda | 0.19.4 |
| foreach | 1.5.2 |
| doParallel | 1.0.17 |
| doRNG | 1.8.3 |
| subplex | 1.8 |
| tidyverse | 1.3.2 |
Arulampalam MS, Maskell S, Gordon N, Clapp T (2002). “A Tutorial on Particle Filters for Online Nonlinear, Non-Gaussian Bayesian Tracking.” IEEE Trans Signal Process, 50, 174–188. https://doi.org/10.1109/78.978374.
Ionides EL, Bretó C, King AA (2006). “Inference for Nonlinear Dynamical Systems.” Proc Natl Acad Sci, 103(49), 18438–18443. https://doi.org/10.1073/pnas.0603181103.
Ionides EL, Nguyen D, Atchadé Y, Stoev S, King AA (2015). “Inference for Dynamic and Latent Variable Models via Iterated, Perturbed Bayes Maps.” Proc Natl Acad Sci, 112(3), 719–724. https://doi.org/10.1073/pnas.1410597112.
McCleery RH, Perrins CM (1991). “Effects of Predation on the Numbers of Great Tits Parus Major.” In Bird population studies. Relevence to conservation and management 129–147. Oxford University Press, Oxford.
Ricker WE (1954). “Stock and Recruitment.” Journal of the Fisheries Research Board of Canada, 11, 559–623. https://doi.org/10.1139/f54-039.